The smuthi.postprocessing package

postprocessing

postprocessing.far_field

Manage post processing steps to evaluate the scattered far field

class smuthi.postprocessing.far_field.FarField(polar_angles='default', azimuthal_angles='default', angular_resolution=None, signal_type='intensity', reference_point=None)

Represent the far field amplitude and far field intensity of an electromagnetic field.

The electric field amplitude \(\mathbf{A}(\theta,\phi)\) is defined by

\[\mathbf{E}(\mathbf{r}) \approx \frac{e^{ikr}}{-ikr} \mathbf{A}(\theta,\phi)\]

for \(kr\rightarrow\infty\), compare equation (3.10) of Bohren and Huffman’s textbook on light scattering.

In the above, \(\mathbf{A}(\theta,\phi)\) is a complex, vector valued function of polar and azimuthal angle. It contains information on the amplitude and phase of the scattered electric field in far field domain.

The intensity \(I_{\Omega,j}(\beta, \alpha)\) is defined by

\[P = \sum_{j=1}^2 \iint \mathrm{d}^2 \Omega \, I_{\Omega,j}(\beta, \alpha),\]

where \(P\) is the radiative power, \(j\) indicates the polarization and \(\mathrm{d}^2 \Omega = \mathrm{d}\alpha \sin\beta \mathrm{d}\beta\) denotes the infinitesimal solid angle.

Parameters:
  • polar_angles (numpy.ndarray) – array of polar angles for plane wave expansions. If ‘default’, use smuthi.fields.default_polar_angles
  • azimuthal_angles (ndarray or str) – array of azimuthal angles for plane wave expansions. If ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • signal_type (str) – Use this field to describe the physical meaning of the power related signal (e.g., ‘intensity’ for standard power flux far fields).
  • reference_point (list or tuple) – [x, y, z]-coordinates of point relative to which the far field is defined
alpha_grid()
Returns:Meshgrid with \(\alpha\) values.
append(other)

Combine two FarField objects with disjoint angular ranges. The other far field is appended to this one.

Parameters:other (FarField) – far field to append to this one.
azimuthal_integral()

Far field intensity as a function of the polar angle cosine only.

\[P = \sum_{j=1}^2 \int \mathrm{d} \cos\beta \, I_{\cos\beta,j}(\beta),\]

with

\[I_{\beta,j}(\beta) = \int \mathrm{d} \alpha \, I_j(\beta, \alpha),\]
Returns:\(I_{\cos\beta,j}(\beta)\) as numpy ndarray. First index is polarization, second is polar angle.
azimuthal_integral_times_sin_beta()

Far field intensity as a function of polar angle only.

\[P = \sum_{j=1}^2 \int \mathrm{d} \beta \, I_{\beta,j}(\beta),\]

with

\[I_{\beta,j}(\beta) = \int \mathrm{d} \alpha \, \sin\beta I_j(\beta, \alpha),\]
Returns:\(I_{\beta,j}(\beta)\) as numpy ndarray. First index is polarization, second is polar angle.
beta_grid()
Returns:Meshgrid with \(\beta\) values.
bottom()

Split far field into top and bottom part.

Returns:FarField object with only the intensity for bottom hemisphere (\(\beta\geq\pi/2\))
electric_field_amplitude()

Evaluate electric field amplitude vector

Returns:Tuple of (A_x, A_y, A_z) numpy.ndarray objects with the Cartesian coordinates of complex electric field amplitude.
integral()

Integrate intensity to obtain total power \(P\).

Returns:\(P_j\) as numpy 1D-array with length 2, the index referring to polarization.
top()

Split far field into top and bottom part.

Returns:FarField object with only the intensity for top hemisphere (\(\beta\leq\pi/2\))
smuthi.postprocessing.far_field.extinction_cross_section(simulation=None, initial_field=None, particle_list=None, layer_system=None, only_l=None, only_m=None, only_pol=None, only_tau=None, extinction_direction='both')

Evaluate the extinction cross section.

Parameters:
  • simulation (smuthi.Simulation.simulation) – Simulation object (optional)
  • initial_field (smuthi.initial_field.PlaneWave) – Plane wave object (optional)
  • particle_list (list) – List of smuthi.particles.Particle objects (optional)
  • layer_system (smuthi.layers.LayerSystem) – Representing the stratified medium
  • only_pol (int) – if set to 0 or 1, only this plane wave polarization (0=TE, 1=TM) is considered
  • only_tau (int) – if set to 0 or 1, only this spherical vector wave polarization (0 — magnetic, 1 — electric) is considered
  • only_l (int) – if set to positive number, only this multipole degree is considered
  • only_m (int) – if set non-negative number, only this multipole order is considered
  • extinction_direction (string) – if set to ‘both’: return full excinction, if to ‘reflection’: extinction of reflected wave, if to ‘transmission’: extinction of transmitted wave. See section on Extinction cross section for details.
Returns:

Extinction cross section.

smuthi.postprocessing.far_field.pwe_to_ff_conversion(vacuum_wavelength, plane_wave_expansion)

Compute the far field of a plane wave expansion object.

Parameters:
  • vacuum_wavelength (float) – Vacuum wavelength in length units.
  • plane_wave_expansion (PlaneWaveExpansion) – Plane wave expansion to convert into far field object.
Returns:

A FarField object containing the far field intensity.

smuthi.postprocessing.far_field.scattered_far_field(vacuum_wavelength, particle_list, layer_system, polar_angles='default', azimuthal_angles='default', angular_resolution=None, reference_point=None)

Evaluate the scattered far field.

Parameters:
  • vacuum_wavelength (float) – in length units
  • particle_list (list) – list of smuthi.Particle objects
  • layer_system (smuthi.layers.LayerSystem) – represents the stratified medium
  • polar_angles (numpy.ndarray or str) – polar angles values (radian). if ‘default’, use smuthi.fields.default_polar_angles
  • azimuthal_angles (numpy.ndarray or str) – azimuthal angle values (radian) if ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • reference_point (list or tuple) – If set to a value other than None, the far field will be calculated with this as its reference point.
Returns:

A smuthi.field_expansion.FarField object of the scattered field.

smuthi.postprocessing.far_field.scattering_cross_section(initial_field, particle_list, layer_system, polar_angles='default', azimuthal_angles='default', angular_resolution=None)

Evaluate and display the differential scattering cross section as a function of solid angle.

Parameters:
  • initial_field (smuthi.initial.PlaneWave) – Initial Plane wave
  • particle_list (list) – scattering particles
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
  • polar_angles (numpy.ndarray or str) – polar angles values (radian). if ‘default’, use smuthi.fields.default_polar_angles
  • azimuthal_angles (numpy.ndarray or str) – azimuthal angle values (radian) if ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
Returns:

A smuthi.field_expansion.FarField object.

smuthi.postprocessing.far_field.total_far_field(initial_field, particle_list, layer_system, polar_angles='default', azimuthal_angles='default', angular_resolution=None, reference_point=None)

Evaluate the total far field, the initial far field and the scattered far field. Cannot be used if initial field is a plane wave.

Parameters:
  • initial_field (smuthi.initial_field.InitialField) – represents the initial field
  • particle_list (list) – list of smuthi.Particle objects
  • layer_system (smuthi.layers.LayerSystem) – represents the stratified medium
  • polar_angles (numpy.ndarray or str) – polar angles values (radian). if ‘default’, use smuthi.fields.default_polar_angles
  • azimuthal_angles (numpy.ndarray or str) – azimuthal angle values (radian) if ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • reference_point (list or tuple) – If set to a value other than None, the far field will be calculated with this as its reference point.
Returns:

A tuple of three smuthi.field_expansion.FarField objects for total, initial and scattered far field. Mind that the scattered far field has no physical meaning and is for illustration purposes only.

smuthi.postprocessing.far_field.total_scattering_cross_section(simulation=None, initial_field=None, particle_list=None, layer_system=None, polar_angles='default', azimuthal_angles='default', angular_resolution=None)

Evaluate the total scattering cross section.

Parameters:
  • simulation (smuthi.Simulation.simulation) – Simulation object (optional)
  • initial_field (smuthi.initial_field.PlaneWave) – Initial Plane wave (optional)
  • particle_list (list) – scattering particles (optional)
  • layer_system (smuthi.layers.LayerSystem) – stratified medium (optional)
  • polar_angles (numpy.ndarray or str) – polar angles values (radian, default None). If None, use smuthi.fields.default_polar_angles
  • azimuthal_angles (numpy.ndarray or str) – azimuthal angle values (radian, default None). If None, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
Returns:

A tuple of smuthi.field_expansion.FarField objects, one for forward scattering (i.e., into the top hemisphere) and one for backward scattering (bottom hemisphere).

postprocessing.graphical_output

Functions to generate plots and animations.

smuthi.postprocessing.graphical_output.compute_near_field(simulation=None, X=None, Y=None, Z=None, type='scatt', chunksize=None, k_parallel='default', azimuthal_angles='default', angular_resolution=None)

Compute a certain component of the electric near field

smuthi.postprocessing.graphical_output.plot_layer_interfaces(dim1min, dim1max, layer_system)

Add lines to plot to display layer system interfaces

Parameters:
  • dim1min (float) – From what x-value plot line
  • dim1max (float) – To what x-value plot line
  • layer_system (smuthi.layers.LayerSystem) – Stratified medium
smuthi.postprocessing.graphical_output.plot_particles(xmin, xmax, ymin, ymax, zmin, zmax, particle_list, draw_circumscribing_sphere, fill_particle=True)

Add circles, ellipses and rectangles to plot to display spheres, spheroids and cylinders.

Parameters:
  • xmin (float) – Minimal x-value of plot
  • xmax (float) – Maximal x-value of plot
  • ymin (float) – Minimal y-value of plot
  • ymax (float) – Maximal y-value of plot
  • zmin (float) – Minimal z-value of plot
  • zmax (float) – Maximal z-value of plot
  • particle_list (list) – List of smuthi.particles.Particle objects
  • draw_circumscribing_sphere (bool) – If true (default), draw a circle indicating the circumscribing sphere of particles.
  • fill_particle (bool) – If true, draw opaque particles.
smuthi.postprocessing.graphical_output.show_far_field(far_field, show_plots=True, show_opts=[{'label': 'far_field'}], save_plots=False, save_opts=None, save_data=False, data_format='hdf5', outputdir='.', flip_downward=True, split=True, log_scale=False)

Display and export the far field.

Parameters:
  • far_field (smuthi.field_expansion.FarField) – Far-field object to show and export
  • show_plots (bool) – Display plots if True
  • show_opts (dict list) – List of dictionaries containing options to be passed to pcolormesh for plotting. If save_plots=True, a 1:1 correspondence between show_opts and save_opts dictionaries is assumed. For simplicity, one can also provide a single show_opts entry that will be applied to all save_opts. The following keys are available (see matplotlib.pyplot.pcolormesh documentation): ‘alpha’ (None) ‘cmap’ (‘inferno’) ‘norm’ (None), is set to matplotlib.colors.LogNorm() if log_scale is True ‘vmin’ (None), applies only to 2D plots ‘vmax’ (None), applies only to 2D plots ‘shading’ (‘nearest’), applies only to 2D plots. ‘gouraud’ is also available ‘linewidth’ (None), applies only to 1D plots ‘linestyle’ (None), applies only to 1D plots ‘marker’ (None), applies only to 1D plots An optional extra key called ‘label’ of type string is shown in the plot title and appended to the associated file if save_plots is True Finally, an optional ‘figsize’ key is available to set the width and height of the figure window (see matplotlib.pyplot.figure documentation)
  • save_plots (bool) – If True, plots are exported to file.
  • save_opts (dict list) – List of dictionaries containing options to be passed to savefig. A 1:1 correspondence between save_opts and show_opts dictionaries is assumed. For simplicity, one can also provide a single save_opts entry that will be applied to all show_opts. The following keys are made available (see matplotlib.pyplot.savefig documentation): ‘dpi’ (None) ‘orientation’ (None) ‘format’ (‘png’), also available: eps, jpeg, jpg, pdf, ps, svg, tif, tiff … ‘transparent’ (False) ‘bbox_inches’ (‘tight’) ‘pad_inches’ (0.1)
  • save_data (bool) – If True, raw data are exported to file
  • data_format (str) – Output data format string, ‘hdf5’ and ‘ascii’ formats are available
  • outputdir (str) – Path to the directory where files are to be saved
  • flip_downward (bool) – If True, represent downward directions as 0-90 deg instead of 90-180
  • split (bool) – If True, show two different plots for upward and downward directions
  • log_scale (bool) – If True, set a logarithmic scale
smuthi.postprocessing.graphical_output.show_near_field(simulation=None, quantities_to_plot=None, show_plots=True, show_opts=None, save_plots=False, save_opts=None, save_data=False, data_format='hdf5', outputdir='.', xmin=0, xmax=0, ymin=0, ymax=0, zmin=0, zmax=0, resolution_step=25, k_parallel='default', azimuthal_angles='default', angular_resolution=None, draw_circumscribing_sphere=True, show_internal_field=False)

Plot the electric near field along a plane. To plot along the xy-plane, specify zmin=zmax and so on.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • quantities_to_plot

    List of strings that specify what to plot. Select from ‘E_x’, ‘E_y’, ‘E_z’, ‘norm(E)’ The list may contain one or more of the following strings:

    ’E_x’ real part of x-component of complex total electric field ‘E_y’ real part of y-component of complex total electric field ‘E_z’ real part of z-component of complex total electric field ‘norm(E)’ norm of complex total electric field

    ’E_scat_x’ real part of x-component of complex scattered electric field ‘E_scat_y’ real part of y-component of complex scattered electric field ‘E_scat_z’ real part of z-component of complex scattered electric field ‘norm(E_scat)’ norm of complex scattered electric field

    ’E_init_x’ real part of x-component of complex initial electric field ‘E_init_y’ real part of y-component of complex initial electric field ‘E_init_z’ real part of z-component of complex initial electric field ‘norm(E_init)’ norm of complex initial electric field

  • show_plots (logical) – If True, plots are shown
  • show_opts (dict list) – List of dictionaries containing options to be passed to imshow for plotting. For each entry in quantities_to_plot, all show_opts dictionaries will be applied. If save_plots=True, a 1:1 correspondence between show_opts and save_opts dictionaries is assumed. For simplicity, one can also provide a single show_opts entry that will be applied to all save_opts. The following keys are made available (see matplotlib.pyplot.imshow documentation): ‘cmap’ defaults to ‘inferno’ for norm quantities and ‘RdYlBu’ otherwise ‘norm’ (None). If a norm is provided, its vmin and vmax take precedence ‘aspect’ (‘equal’) ‘interpolation’ (None), also available: bilinear, bicubic, spline16, quadric, … ‘alpha’ (None) ‘vmin’ (None), will be set to 0 for norm quantities and -vmax otherwise ‘vmax’ initialized with the max of the quantity to plot ‘origin’ (‘lower’) ‘extent’ calculated automatically based on plotting coordinate limits An optional extra key called ‘label’ of type string is shown in the plot title and appended to the associated file if save_plots is True Finally, an optional ‘figsize’ key is available to set the width and height of the figure window (see matplotlib.pyplot.figure documentation)
  • save_plots (logical) – If True, plots are exported to file.
  • save_opts (dict list) – List of dictionaries containing options to be passed to savefig. For each entry in quantities_to_plot, all save_opts dictionaries will be applied. A 1:1 correspondence between save_opts and show_opts dictionaries is assumed. For simplicity, one can also provide a single save_opts entry that will be applied to all show_opts. The following keys are made available (see matplotlib.pyplot.savefig documentation): ‘dpi’ (None) ‘orientation’ (None) ‘format’ (‘png’), also available: eps, jpeg, jpg, pdf, ps, svg, tif, tiff … ‘transparent’ (False) ‘bbox_inches’ (‘tight’) ‘pad_inches’ (0.1) Passing ‘gif’ as one of the format values will result in an animation if the quantity to plot is of non-norm type
  • save_data (logical) – If True, raw data are exported to file
  • data_format (str) – Output data format string, ‘hdf5’ and ‘ascii’ formats are available
  • outputdir (str) – Path to directory where to save the export files
  • xmin (float) – Plot from that x (length unit)
  • xmax (float) – Plot up to that x (length unit)
  • ymin (float) – Plot from that y (length unit)
  • ymax (float) – Plot up to that y (length unit)
  • zmin (float) – Plot from that z (length unit)
  • zmax (float) – Plot up to that z (length unit)
  • resolution_step (float) – Compute the field with that spatial resolution (length unit, distance between computed points), can be a tuple for [resx, resy, resz]
  • k_parallel (numpy.ndarray or str) – in-plane wavenumbers for the plane wave expansion if ‘default’, use smuthi.fields.default_Sommerfeld_k_parallel_array
  • azimuthal_angles (numpy.ndarray or str) – azimuthal angles for the plane wave expansion if ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • draw_circumscribing_sphere (bool) – If true (default), draw a circle indicating the circumscribing sphere of particles.
  • show_internal_field (bool) – If true, compute also the field inside the particles (only for spheres)
smuthi.postprocessing.graphical_output.show_scattered_far_field(simulation, show_plots=True, show_opts=[{'label': 'scattered_far_field'}], save_plots=False, save_opts=None, save_data=False, data_format='hdf5', outputdir='.', flip_downward=True, split=True, log_scale=False, polar_angles='default', azimuthal_angles='default', angular_resolution=None)

Display and export the scattered far field.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • show_plots (bool) – Display plots if True
  • show_opts (dict list) – List of dictionaries containing options to be passed to pcolormesh for plotting. If save_plots=True, a 1:1 correspondence between show_opts and save_opts dictionaries is assumed. For simplicity, one can also provide a single show_opts entry that will be applied to all save_opts. The following keys are available (see matplotlib.pyplot.pcolormesh documentation): ‘alpha’ (None) ‘cmap’ (‘inferno’) ‘norm’ (None), is set to matplotlib.colors.LogNorm() if log_scale is True ‘vmin’ (None), applies only to 2D plots ‘vmax’ (None), applies only to 2D plots ‘shading’ (‘nearest’), applies only to 2D plots. ‘gouraud’ is also available ‘linewidth’ (None), applies only to 1D plots ‘linestyle’ (None), applies only to 1D plots ‘marker’ (None), applies only to 1D plots An optional extra key called ‘label’ of type string is shown in the plot title and appended to the associated file if save_plots is True
  • save_plots (bool) – If True, plots are exported to file.
  • save_opts (dict list) – List of dictionaries containing options to be passed to savefig. A 1:1 correspondence between save_opts and show_opts dictionaries is assumed. For simplicity, one can also provide a single save_opts entry that will be applied to all show_opts. The following keys are made available (see matplotlib.pyplot.savefig documentation): ‘dpi’ (None) ‘orientation’ (None) ‘format’ (‘png’), also available: eps, jpeg, jpg, pdf, ps, svg, tif, tiff … ‘transparent’ (False) ‘bbox_inches’ (‘tight’) ‘pad_inches’ (0.1)
  • save_data (bool) – If True, raw data are exported to file
  • data_format (str) – Output data format string, ‘hdf5’ and ‘ascii’ formats are available
  • outputdir (str) – Path to the directory where files are to be saved
  • flip_downward (bool) – If True, represent downward directions as 0-90 deg instead of 90-180
  • split (bool) – If True, show two different plots for upward and downward directions
  • log_scale (bool) – If True, set a logarithmic scale
  • polar_angles (numpy.ndarray or str) – Polar angles values (radian). If ‘default’, use smuthi.fields.default_polar_angles
  • azimuthal_angles (numpy.ndarray or str) – Azimuthal angle values (radian). If ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
smuthi.postprocessing.graphical_output.show_scattering_cross_section(simulation, show_plots=True, show_opts=[{'label': 'scattering_cross_section'}], save_plots=False, save_opts=None, save_data=False, data_format='hdf5', outputdir='.', flip_downward=True, split=True, log_scale=False, polar_angles='default', azimuthal_angles='default', angular_resolution=None)

Display and export the differential scattering cross section.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • show_plots (bool) – Display plots if True
  • show_opts (dict list) – List of dictionaries containing options to be passed to pcolormesh for plotting. If save_plots=True, a 1:1 correspondence between show_opts and save_opts dictionaries is assumed. For simplicity, one can also provide a single show_opts entry that will be applied to all save_opts. The following keys are available (see matplotlib.pyplot.pcolormesh documentation): ‘alpha’ (None) ‘cmap’ (‘inferno’) ‘norm’ (None), is set to matplotlib.colors.LogNorm() if log_scale is True ‘vmin’ (None), applies only to 2D plots ‘vmax’ (None), applies only to 2D plots ‘shading’ (‘nearest’), applies only to 2D plots. ‘gouraud’ is also available ‘linewidth’ (None), applies only to 1D plots ‘linestyle’ (None), applies only to 1D plots ‘marker’ (None), applies only to 1D plots An optional extra key called ‘label’ of type string is shown in the plot title and appended to the associated file if save_plots is True
  • save_plots (bool) – If True, plots are exported to file.
  • save_opts (dict list) – List of dictionaries containing options to be passed to savefig. A 1:1 correspondence between save_opts and show_opts dictionaries is assumed. For simplicity, one can also provide a single save_opts entry that will be applied to all show_opts. The following keys are made available (see matplotlib.pyplot.savefig documentation): ‘dpi’ (None) ‘orientation’ (None) ‘format’ (‘png’), also available: eps, jpeg, jpg, pdf, ps, svg, tif, tiff … ‘transparent’ (False) ‘bbox_inches’ (‘tight’) ‘pad_inches’ (0.1)
  • save_data (bool) – If True, raw data are exported to file
  • data_format (str) – Output data format string, ‘hdf5’ and ‘ascii’ formats are available
  • outputdir (str) – Path to the directory where files are to be saved
  • flip_downward (bool) – If True, represent downward directions as 0-90 deg instead of 90-180
  • split (bool) – If True, show two different plots for upward and downward directions
  • log_scale (bool) – If True, set a logarithmic scale
  • polar_angles (numpy.ndarray or str) – Polar angles values (radian). If ‘default’, use smuthi.fields.default_polar_angles
  • azimuthal_angles (numpy.ndarray or str) – Azimuthal angle values (radian). If ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
smuthi.postprocessing.graphical_output.show_total_far_field(simulation, show_plots=True, show_opts=[{'label': 'total_far_field'}], save_plots=False, save_opts=None, save_data=False, data_format='hdf5', outputdir='.', flip_downward=True, split=True, log_scale=False, polar_angles='default', azimuthal_angles='default', angular_resolution=None)

Display and export the total far field. This function cannot be used if the inital field is a plane wave.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • show_plots (bool) – Display plots if True
  • show_opts (dict list) – List of dictionaries containing options to be passed to pcolormesh for plotting. If save_plots=True, a 1:1 correspondence between show_opts and save_opts dictionaries is assumed. For simplicity, one can also provide a single show_opts entry that will be applied to all save_opts. The following keys are available (see matplotlib.pyplot.pcolormesh documentation): ‘alpha’ (None) ‘cmap’ (‘inferno’) ‘norm’ (None), is set to matplotlib.colors.LogNorm() if log_scale is True ‘vmin’ (None), applies only to 2D plots ‘vmax’ (None), applies only to 2D plots ‘shading’ (‘nearest’), applies only to 2D plots. ‘gouraud’ is also available ‘linewidth’ (None), applies only to 1D plots ‘linestyle’ (None), applies only to 1D plots ‘marker’ (None), applies only to 1D plots An optional extra key called ‘label’ of type string is shown in the plot title and appended to the associated file if save_plots is True
  • save_plots (bool) – If True, plots are exported to file.
  • save_opts (dict list) – List of dictionaries containing options to be passed to savefig. A 1:1 correspondence between save_opts and show_opts dictionaries is assumed. For simplicity, one can also provide a single save_opts entry that will be applied to all show_opts. The following keys are made available (see matplotlib.pyplot.savefig documentation): ‘dpi’ (None) ‘orientation’ (None) ‘format’ (‘png’), also available: eps, jpeg, jpg, pdf, ps, svg, tif, tiff … ‘transparent’ (False) ‘bbox_inches’ (‘tight’) ‘pad_inches’ (0.1)
  • save_data (bool) – If True, raw data are exported to file
  • data_format (str) – Output data format string, ‘hdf5’ and ‘ascii’ formats are available
  • outputdir (str) – Path to the directory where files are to be saved
  • flip_downward (bool) – If True, represent downward directions as 0-90 deg instead of 90-180
  • split (bool) – If True, show two different plots for upward and downward directions
  • log_scale (bool) – If True, set a logarithmic scale
  • polar_angles (numpy.ndarray or str) – Polar angles values (radian). If ‘default’, use smuthi.fields.default_polar_angles
  • azimuthal_angles (numpy.ndarray or str) – Azimuthal angle values (radian). If ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range

postprocessing.internal_field


Manage post processing steps to evaluate the electric field inside a sphere

smuthi.postprocessing.internal_field.internal_field_piecewise_expansion(vacuum_wavelength, particle_list, layer_system)

Compute a piecewise field expansion of the internal field of spheres.

Parameters:
  • vacuum_wavelength (float) – vacuum wavelength
  • particle_list (list) – list of smuthi.particles.Particle objects
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
Returns:

internal field as smuthi.field_expansion.PiecewiseFieldExpansion object

postprocessing.scattered_field

Manage post processing steps to evaluate the scattered electric field

smuthi.postprocessing.scattered_field.evaluate_scattered_field_stat_phase_approx(x, y, z, vacuum_wavelength, particle_list, layer_system)

Evaluate the scattered electric field for N particles on a substrate. The substrate reflection is evaluated by means of the stationary phase approximation, as presented in “A quick way to approximate a Sommerfeld-Weyl_type Sommerfeld integral” by W.C. Chew (1988).

See also the technical note “Usage of the stationary phase approximation in SMUTHI” by A. Egel (2020)

The stationary phase approximation is expected to yield good results for field points far away from the particles.

Note: This function assumes that the particles are located in the upper layer of a two-layer system (particles on substrate). For other cases, this function does not apply. ****************************************************************************************************************

Parameters:
  • x (float or numpy.ndarray) – x-coordinates of query points
  • y (float or numpy.ndarray) – y-coordinates of query points
  • z (float or numpy.ndarray) – z-coordinates of query points
  • vacuum_wavelength (float) – Vacuum wavelength \(\lambda\) (length unit)
  • particle_list (list) – List of Particle objects
  • layer_system (smuthi.layers.LayerSystem) – Stratified medium
Returns:

Tuple of (E_x, E_y, E_z) numpy.ndarray objects with the Cartesian coordinates of complex electric field.

smuthi.postprocessing.scattered_field.scattered_field_piecewise_expansion(vacuum_wavelength, particle_list, layer_system, k_parallel='default', azimuthal_angles='default', angular_resolution=None, layer_numbers=None)

Compute a piecewise field expansion of the scattered field.

Parameters:
  • vacuum_wavelength (float) – vacuum wavelength
  • particle_list (list) – list of smuthi.particles.Particle objects
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
  • k_parallel (numpy.ndarray or str) – in-plane wavenumbers array. if ‘default’, use smuthi.fields.default_Sommerfeld_k_parallel_array
  • azimuthal_angles (numpy.ndarray or str) – azimuthal angles array if ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • layer_numbers (list) – if specified, append only plane wave expansions for these layers
Returns:

scattered field as smuthi.field_expansion.PiecewiseFieldExpansion object

smuthi.postprocessing.scattered_field.scattered_field_pwe(vacuum_wavelength, particle_list, layer_system, layer_number, k_parallel='default', azimuthal_angles='default', angular_resolution=None, include_direct=True, include_layer_response=True, only_l=None, only_m=None, only_pol=None, only_tau=None)

Calculate the plane wave expansion of the scattered field of a set of particles.

Parameters:
  • vacuum_wavelength (float) – Vacuum wavelength (length unit)
  • particle_list (list) – List of Particle objects
  • layer_system (smuthi.layers.LayerSystem) – Stratified medium
  • layer_number (int) – Layer number in which the plane wave expansion should be valid
  • k_parallel (numpy.ndarray or str) – in-plane wavenumbers array. if ‘default’, use smuthi.fields.default_Sommerfeld_k_parallel_array
  • azimuthal_angles (numpy.ndarray or str) – azimuthal angles array if ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • include_direct (bool) – If True, include the direct scattered field
  • include_layer_response (bool) – If True, include the layer system response
  • only_pol (int) – if set to 0 or 1, only this plane wave polarization (0=TE, 1=TM) is considered
  • only_tau (int) – if set to 0 or 1, only this spherical vector wave polarization (0 — magnetic, 1 — electric) is considered
  • only_l (int) – if set to positive number, only this multipole degree is considered
  • only_m (int) – if set to non-negative number, only this multipole order is considered
Returns:

A tuple of PlaneWaveExpansion objects for upgoing and downgoing waves.

postprocessing.power_flux

Manage post processing steps to evaluate power flux

smuthi.postprocessing.power_flux.power_flux_through_zplane(vacuum_wavelength, z, upgoing_pwe=None, downgoing_pwe=None)

Evaluate time averaged power flux though a plane of z=const.

Parameters:
  • vacuum_wavelength (float) – Vacuum wavelength in length units.
  • z (float) – plane height z
  • upgoing_pwe (PlaneWaveExpansion) – of kind “upgoing”
  • downgoing_pwe (PlaneWaveExpansion) – of kind “downgoing”
Returns:

Time averaged energy flux.